dspawpy.diffusion package¶
Submodules¶
dspawpy.diffusion.neb module¶
- class dspawpy.diffusion.neb.NEB(initial_structure, final_structure, nimages)¶
Bases:
object
NEB插值算法
- Parameters:
initial_structure (Structure) -- 初态
final_structure (Structure) -- 终态
nimages (int) -- 中间构型数
- idpp_interpolate(maxiter=1000, tol=1e-05, gtol=0.001, step_size=0.05, max_disp=0.05, spring_const=5.0)¶
- linear_interpolate()¶
- dspawpy.diffusion.neb.write_neb_structures(structures, coords_are_cartesian=True, fmt: str = 'as', path: str = '.', prefix='structure')¶
插值并生成中间构型文件
- Parameters:
structures (list) -- 构型列表
coords_are_cartesian (bool) -- 坐标是否为笛卡尔坐标
fmt (str) -- 结构文件类型,默认为 "as"
path (str) -- 保存路径
prefix (str) -- 文件名前缀,默认为 "structure",这样的话生成的就是 structure00.as, structure01.as, ...
- Returns:
保存构型文件
- Return type:
file
Examples
先读取as文件创建structure对象
>>> from dspawpy.io.structure import read >>> init_struct = read("/data/home/hzw1002/dspawpy_repo/test/2.15/00/structure00.as")[0] >>> final_struct = read("/data/home/hzw1002/dspawpy_repo/test/2.15/04/structure04.as")[0]
然后,插值并生成中间构型文件
>>> from dspawpy.diffusion.neb import NEB,write_neb_structures >>> neb = NEB(init_struct,final_struct,8) >>> structures = neb.linear_interpolate() #线性插值
插值完成的构型可指定保存到neb文件夹下
>>> write_neb_structures(structures, fmt="as", path="/data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures") ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/00/structure00.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/01/structure01.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/02/structure02.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/03/structure03.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/04/structure04.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/05/structure05.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/06/structure06.as ==> /data/home/hzw1002/dspawpy_repo/test/out/11neb_interpolate_structures/07/structure07.as
dspawpy.diffusion.nebtools module¶
- dspawpy.diffusion.nebtools.get_distance(spo1: ndarray, spo2: ndarray, lat1: ndarray, lat2: ndarray)¶
根据两个结构的分数坐标和晶胞计算距离
- Parameters:
spo1 (np.ndarray) -- 分数坐标列表1
spo2 (np.ndarray) -- 分数坐标列表2
lat1 (np.ndarray) -- 晶胞1
lat2 (np.ndarray) -- 晶胞2
- Returns:
距离
- Return type:
float
Examples
先读取结构信息
>>> from dspawpy.io.structure import read >>> s1 = read('/data/home/hzw1002/dspawpy_repo/test/2.15/01/structure01.as')[0] >>> s2 = read('/data/home/hzw1002/dspawpy_repo/test/2.15/02/structure02.as')[0]
计算两个构型的距离
>>> from dspawpy.diffusion.nebtools import get_distance >>> dist = get_distance(s1.frac_coords, s2.frac_coords, s1.lattice.matrix, s2.lattice.matrix) >>> print('两个构型的距离为:', dist, 'Angstrom') 两个构型的距离为: 0.476972808803491 Angstrom
- dspawpy.diffusion.nebtools.get_neb_subfolders(directory: str = '.', return_abs=False)¶
将directory路径下的子文件夹名称列表按照数字大小排序
仅保留形如00,01数字类型的NEB子文件夹路径
- Parameters:
subfolders (list) -- 子文件夹名称列表
return_abs (bool, optional) -- 是否返回绝对路径, 默认否
- Returns:
subfolders -- 排序后的子文件夹名称列表
- Return type:
list
Examples
>>> from dspawpy.diffusion.nebtools import get_neb_subfolders >>> directory = '/data/home/hzw1002/dspawpy_repo/test/2.15' >>> get_neb_subfolders(directory) ['00', '01', '02', '03', '04']
- dspawpy.diffusion.nebtools.plot_barrier(datafile: str = 'neb.h5', directory: str = None, ri: float = None, rf: float = None, ei: float = None, ef: float = None, method: str = 'PchipInterpolator', figname: str = 'neb_barrier.png', show: bool = True, raw=False, **kwargs)¶
调用 scipy.interpolate 插值算法,拟合NEB能垒并绘图
- Parameters:
datafile (str) -- neb.h5或neb.json文件路径
directory (str) -- NEB计算路径
ri (float) -- 初态反应坐标
rf (float) -- 末态反应坐标
ei (float) -- 初态自洽能量
ef (float) -- 末态自洽能量
method (str, optional) -- 插值算法, 默认'PchipInterpolator'
figname (str, optional) -- 能垒图名称, 默认'neb_barrier.png'
show (bool, optional) -- 是否展示交互界面, 默认True
raw (bool, optional) -- 是否返回绘图数据到csv
- Raises:
ImportError -- 指定了scipy.interpolate中不存在的插值算法
ValueError -- 传递给插值算法的参数不符合该算法要求
Examples
>>> from dspawpy.diffusion.nebtools import plot_barrier >>> import matplotlib.pyplot as plt
对比不同插值算法
>>> plot_barrier(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', method='interp1d', kind=2, figname=None, show=False) >>> plot_barrier(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', method='interp1d', kind=3, figname=None, show=False) >>> plot_barrier(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', method='CubicSpline', figname=None, show=False) >>> plot_barrier(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', method='pchip', figname='../APIout/barrier_comparison.png', show=False) ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/barrier_comparison.png
尝试读取neb.h5文件或neb.json文件
>>> plot_barrier(datafile='/data/home/hzw1002/dspawpy_repo/test/2.15/neb.h5', method='pchip', figname='../APIout/barrier_h5.png', show=False) ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/barrier_h5.png >>> plot_barrier(datafile='/data/home/hzw1002/dspawpy_repo/test/2.15/neb.json', method='pchip', figname='../APIout/barrier_json.png', show=False) ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/barrier_json.png
- dspawpy.diffusion.nebtools.plot_barrier_old(datafile: str = 'neb.h5', directory: str = None, ri: float = None, rf: float = None, ei: float = None, ef: float = None, method: str = 'PchipInterpolator', figname: str = 'neb_barrier.png', show: bool = True, raw=False, **kwargs)¶
与plot_barrier()唯一的区别在于反应坐标视为累加值,而不是单个值
- dspawpy.diffusion.nebtools.plot_neb_converge(neb_dir: str, image_key: str = '01', show: bool = True, figname: str = 'neb_conv.png', raw=False)¶
指定NEB计算路径,绘制NEB收敛过程图
- Parameters:
neb_dir (str) -- neb.h5 / neb.json 文件路径或者包含 neb.h5 / neb.json 文件的文件夹路径
image_key (str) -- 第几个构型,默认 "01"
show (bool) -- 是否交互绘图
image_name (str) -- NEB收敛图名称,默认 "neb_conv.png"
raw (bool) -- 是否输出绘图数据到csv文件
- Returns:
ax1, ax2 -- 两个子图的Axes对象
- Return type:
matplotlib.axes.Axes
Examples
>>> from dspawpy.diffusion.nebtools import plot_neb_converge >>> plot_neb_converge(neb_dir='/data/home/hzw1002/dspawpy_repo/test/2.15', image_key='01', figname='/data/home/hzw1002/dspawpy_repo/test/out/neb_converge1.png',show=False) ==> /data/home/hzw1002/dspawpy_repo/test/out/neb_converge1.png (<Axes: xlabel='Number of ionic step', ylabel='Force (eV/Å)'>, <Axes: xlabel='Number of ionic step', ylabel='Energy (eV)'>)
- dspawpy.diffusion.nebtools.printef(directory)¶
打印NEB计算时各构型的能量和受力
- Parameters:
directory (str) -- NEB计算的目录,默认为当前目录
- Return type:
打印各构型的能量和受力
Examples
>>> from dspawpy.diffusion.nebtools import printef >>> printef(directory='/data/home/hzw1002/dspawpy_repo/test/2.15') Force(eV/Å) RC(Å) Energy(eV) E-E0(eV) 00 0.180272 0.000000 -39637.098409 0.000000 01 0.026337 0.542789 -39637.018595 0.079814 02 0.024798 1.086800 -39636.880144 0.218265 03 0.234429 1.588367 -39636.998366 0.100043 04 0.014094 2.089212 -39637.089994 0.008414
- dspawpy.diffusion.nebtools.restart(directory: str = '.', output: str = 'bakfile')¶
将旧NEB任务归档压缩,并在原路径下准备续算
- Parameters:
directory (str) -- 旧NEB任务所在路径,默认当前路径
output (str) -- 备份文件夹路径,默认将在当前路径新建一个bakfile文件夹用于备份; 也可以任意指定一个路径,但不能与当前路径相同
Examples
>>> from dspawpy.diffusion.nebtools import restart >>> restart(directory='/data/home/hzw1002/dspawpy_repo/test/neb_temp_back', output='/data/home/hzw1002/dspawpy_repo/test/out/neb_backup') ==> /data/home/hzw1002/dspawpy_repo/test/out/neb_backup
续算准备工作可能需要较长时间才能完成,请耐心等待
- dspawpy.diffusion.nebtools.set_pbc(spo)¶
根据周期性边界条件将分数坐标分量移入 [-0.5, 0.5) 区间
- Parameters:
spo (np.ndarray or list) -- 分数坐标列表
- Returns:
pbc_spo -- 符合周期性边界条件的分数坐标列表
- Return type:
np.ndarray
Examples
>>> from dspawpy.diffusion.nebtools import set_pbc >>> set_pbc([-0.6, 1.2, 2.3]) array([0.4, 0.2, 0.3])
- dspawpy.diffusion.nebtools.summary(directory: str = '.', raw=False, show_converge=False, outdir: str = None, **kwargs)¶
NEB任务完成总结,依次执行以下步骤:
打印各构型受力、反应坐标、能量、与初始构型的能量差
绘制能垒图
绘制并保存结构优化过程的能量和受力收敛过程图
- Parameters:
directory (str) -- NEB路径, 默认当前路径
raw (bool) -- 是否保存绘图数据到csv文件
show_converge (bool) -- 是否展示结构优化过程的能量和受力收敛过程图,默认不展示
outdir (str) -- 收敛过程图保存路径,默认为directory
**kwargs (dict) -- 传递给plot_barrier的参数
Examples
>>> from dspawpy.diffusion.nebtools import summary >>> directory = '/data/home/hzw1002/dspawpy_repo/test/2.15' # NEB计算路径,默认当前路径 >>> summary(directory, show=False, figname='../APIout/neb_barrier.png') Force(eV/Å) RC(Å) Energy(eV) E-E0(eV) 00 0.180272 0.000000 -39637.098409 0.000000 01 0.026337 0.542789 -39637.018595 0.079814 02 0.024798 1.086800 -39636.880144 0.218265 03 0.234429 1.588367 -39636.998366 0.100043 04 0.014094 2.089212 -39637.089994 0.008414 ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_barrier.png ==> /data/home/hzw1002/dspawpy_repo/test/2.15/01/converge.png... ==> /data/home/hzw1002/dspawpy_repo/test/2.15/01/converge.png ==> /data/home/hzw1002/dspawpy_repo/test/2.15/02/converge.png... ==> /data/home/hzw1002/dspawpy_repo/test/2.15/02/converge.png ==> /data/home/hzw1002/dspawpy_repo/test/2.15/03/converge.png... ==> /data/home/hzw1002/dspawpy_repo/test/2.15/03/converge.png
若inifin=false,用户必须将自洽的scf.h5或system.json放到初末态子文件夹中
- dspawpy.diffusion.nebtools.write_json_chain(preview: bool = False, directory: str = '.', step: int = -1, dst: str = None, ignorels=False)¶
NEB计算或者初始插值后,读取信息,保存为 neb_chain*.json 文件
用 Device Studio 打开该文件可以观察结构等信息
- Parameters:
preview (bool) -- 是否预览模式,默认否
directory (str) -- 计算结果所在目录. 默认当前路径
step (int) --
离子步编号. 默认-1,读取整个NEB计算过程信息;此时优先级为latestStructure.as > h5 > json
0表示初插结构(未完成离子步);
1表示第一个离子步,以此类推
dst (str) -- 保存路径,默认为directory
ignorels (bool) -- 当step=-1时是否忽略latestStructure.as文件,默认否
Examples
>>> from dspawpy.diffusion.nebtools import write_json_chain
NEB计算完成后要观察轨迹变化全过程,只需指定NEB计算路径即可
>>> write_json_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', dst='../APIout') structure info collected from latestStructure.as ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_last.json
NEB计算完成后要观察第n离子步结构,请设置step为n,注意step从1开始计数
>>> write_json_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', step=1, dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_1.json
如果您指定的step数超过NEB实际完成的离子步,将会自动修改为最后一步
>>> write_json_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', step=10, dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_10.json
另外,如需预览初插结构,请将preview设置为True,并将directory指定为NEB计算主路径
>>> write_json_chain(preview=True, directory='/data/home/hzw1002/dspawpy_repo/test/2.15', dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_init.json
- dspawpy.diffusion.nebtools.write_movie_json(preview: bool = False, directory: str = '.', step: int = -1, dst: str = None)¶
- dspawpy.diffusion.nebtools.write_xyz(preview: bool = False, directory: str = '.', step: int = -1, dst: str = None)¶
- dspawpy.diffusion.nebtools.write_xyz_chain(preview: bool = False, directory: str = '.', step: int = -1, dst: str = None, ignorels=False)¶
将NEB结构链条写成xyz轨迹文件用于可视化
- Parameters:
preview (bool) -- 是否预览模式,默认否
directory (str) -- 计算结果所在目录. 默认当前路径
step (int) --
离子步编号. 默认-1,读取整个NEB计算过程信息;此时优先级为latestStructure.as > h5 > json
0表示初插结构(未完成离子步);
1表示第一个离子步,以此类推
dst (str) -- 保存路径,默认为directory
ignorels (bool) -- 当step=-1时是否忽略latestStructure.as文件,默认否
Examples
>>> from dspawpy.diffusion.nebtools import write_xyz_chain
NEB计算完成后要观察轨迹变化全过程,只需指定NEB计算路径即可
>>> write_xyz_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', dst='../APIout') structure info collected from latestStructure.as ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_last.xyz
NEB计算完成后要观察第n离子步结构,请设置step为n,注意step从1开始计数
>>> write_xyz_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', step=1, dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_1.xyz
如果您指定的step数超过NEB实际完成的离子步,将会自动修改为最后一步
>>> write_xyz_chain(directory='/data/home/hzw1002/dspawpy_repo/test/2.15', step=10, dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_10.xyz
另外,如需预览初插结构,请将preview设置为True,并将directory指定为NEB计算主路径
>>> write_xyz_chain(preview=True, directory='/data/home/hzw1002/dspawpy_repo/test/2.15', dst='../APIout') ==> /data/home/hzw1002/dspawpy_repo/dspaw-manual-cn/APIout/neb_chain_init.xyz
dspawpy.diffusion.pathfinder module¶
- class dspawpy.diffusion.pathfinder.IDPPSolver(structures)¶
Bases:
object
A solver using image dependent pair potential (IDPP) algo to get an improved initial NEB path. For more details about this algo, please refer to Smidstrup et al., J. Chem. Phys. 140, 214106 (2014).
- classmethod from_endpoints(endpoints, nimages=5, sort_tol=1.0)¶
A class method that starts with end-point structures instead. The initial guess for the IDPP algo is then constructed using linear interpolation.
- Parameters:
endpoints (list of Structure objects) -- The two end-point structures.
nimages (int) -- Number of images between the two end-points.
sort_tol (float) -- Distance tolerance (in Angstrom) used to match the atomic indices between start and end structures. Need to increase the value in some cases.
- static get_unit_vector(vec)¶
- run(maxiter=1000, tol=1e-05, gtol=0.001, step_size=0.05, max_disp=0.05, spring_const=5.0, species=None)¶
Perform iterative minimization of the set of objective functions in an NEB-like manner. In each iteration, the total force matrix for each image is constructed, which comprises both the spring forces and true forces. For more details about the NEB approach, please see the references, e.g. Henkelman et al., J. Chem. Phys. 113, 9901 (2000).
- Parameters:
maxiter (int) -- Maximum number of iterations in the minimization process.
tol (float) -- Tolerance of the change of objective functions between consecutive steps.
gtol (float) -- Tolerance of maximum force component (absolute value).
step_size (float) -- Step size associated with the displacement of the atoms during the minimization process.
max_disp (float) -- Maximum allowed atomic displacement in each iteration.
spring_const (float) -- A virtual spring constant used in the NEB-like relaxation process that yields so-called IDPP path.
species (list of string) -- If provided, only those given species are allowed to move. The atomic positions of other species are obtained via regular linear interpolation approach.
- Returns:
[Structure] Complete IDPP path (including end-point structures)